Iterative Solvers

4. Iterative Solvers#

So far we have used direct solvers to solve the linear system of equations. Although a direct solver can profit from the sparse matrix, it’s arithmetic complexity is sub-optimal (i.e. not linear in the number of degrees of freedom). For large-scale problems iterative solvers are a must.

The conjugate gradient (cg) method is the standard method for symmetric and positive definite matrices. It’s convergence rate depends on a preconditioner, what is a cheap approximative inverse to the matrix.

from ngsolve import *
from ngsolve.webgui import Draw

We generate a 3D geometry and mesh using the OCC constructive solid geometry (CSG) modeler:

from netgen.occ import *
cube = Box((0,0,0),(1,1,1))
cyl = Cylinder((0,0.5,0.5),X, r=0.2, h=1)
cube.faces.name = "outer"
cyl.faces.name = "cyl"
shape = cube-cyl
Draw(shape);

Generate a mesh, and perform uniform mesh refinement:

ngmesh = OCCGeometry(shape).GenerateMesh(maxh=0.1)
for l in range(1):
    ngmesh.Refine()
mesh = Mesh(ngmesh)
mesh.Curve(3)
Draw (mesh);
fes = H1(mesh, order=3, dirichlet="outer", wb_withedges=False)
print ("we have", fes.ndof, "unknowns")
u,v = fes.TnT()

a = BilinearForm(grad(u)*grad(v)*dx)
f = LinearForm(v*dx)

# c = Preconditioner(a, "direct", inverse="sparsecholesky")
c = Preconditioner(a, "local")
# c = Preconditioner(a, "bddc")

gfu = GridFunction(fes)
we have 177430 unknowns

assemble system and setup preconditioner in parallel:

with TaskManager():
    a.Assemble()
    f.Assemble()

solve the system using the preconditioned conjugate gradient method:

from ngsolve.krylovspace import CGSolver

with TaskManager():
    inv = CGSolver(mat=a.mat, pre=c.mat, printrates='\r', maxiter=400)
    gfu.vec.data = inv * f.vec
CG iteration 1, residual = 0.024273223974312075     
CG iteration 2, residual = 0.037204475895342     
CG iteration 3, residual = 0.048296081732792415     
CG iteration 4, residual = 0.04783145202720861     
CG iteration 5, residual = 0.03503614777663431     
CG iteration 6, residual = 0.029378085186036507     
CG iteration 7, residual = 0.031006371846374867     
CG iteration 8, residual = 0.02881189922538991     
CG iteration 9, residual = 0.02699139348804369     
CG iteration 10, residual = 0.025607813233479118     
CG iteration 11, residual = 0.02305922750715943     
CG iteration 12, residual = 0.0201887142527063     
CG iteration 13, residual = 0.01767785669092636     
CG iteration 14, residual = 0.01570390414750895     
CG iteration 15, residual = 0.014508920088754974     
CG iteration 16, residual = 0.013060236654630367     
CG iteration 17, residual = 0.011487172720698012     
CG iteration 18, residual = 0.010359419640115405     
CG iteration 19, residual = 0.00963928071857335     
CG iteration 20, residual = 0.00898378451986901     
CG iteration 21, residual = 0.00825049617120831     
CG iteration 22, residual = 0.007352872683189113     
CG iteration 23, residual = 0.006439789142169032     
CG iteration 24, residual = 0.005596514201958335     
CG iteration 25, residual = 0.0047662198736664165     
CG iteration 26, residual = 0.004022073286663704     
CG iteration 27, residual = 0.0034196513449910805     
CG iteration 28, residual = 0.002940004822205984     
CG iteration 29, residual = 0.002508333270524199     
CG iteration 30, residual = 0.002104649927409852     
CG iteration 31, residual = 0.0017484379251981052     
CG iteration 32, residual = 0.001439742585061398     
CG iteration 33, residual = 0.0011912144262735038     
CG iteration 34, residual = 0.0009893154642094018     
CG iteration 35, residual = 0.0008329679456139166     
CG iteration 36, residual = 0.0007272861126943781     
CG iteration 37, residual = 0.0006352763333349189     
CG iteration 38, residual = 0.0005645223644237872     
CG iteration 39, residual = 0.0004942564003676567     
CG iteration 40, residual = 0.0004278025775336734     
CG iteration 41, residual = 0.00036683687047396385     
CG iteration 42, residual = 0.0003112297948991816     
CG iteration 43, residual = 0.00026555485540291645     
CG iteration 44, residual = 0.00023057965122062098     
CG iteration 45, residual = 0.0002081072052177469     
CG iteration 46, residual = 0.00019208909401800619     
CG iteration 47, residual = 0.0001785379289683779     
CG iteration 48, residual = 0.00016491550344460986     
CG iteration 49, residual = 0.00014983626453891538     
CG iteration 50, residual = 0.0001338648614332047     
CG iteration 51, residual = 0.0001159879659544772     
CG iteration 52, residual = 9.959344694700486e-05     
CG iteration 53, residual = 8.519641137784852e-05     
CG iteration 54, residual = 7.322459830883046e-05     
CG iteration 55, residual = 6.417417797477809e-05     
CG iteration 56, residual = 5.69968450773375e-05     
CG iteration 57, residual = 5.166607912750669e-05     
CG iteration 58, residual = 4.711071719270879e-05     
CG iteration 59, residual = 4.323161539942041e-05     
CG iteration 60, residual = 3.952383801108811e-05     
CG iteration 61, residual = 3.548362305135954e-05     
CG iteration 62, residual = 3.115776481808804e-05     
CG iteration 63, residual = 2.6858896178966248e-05     
CG iteration 64, residual = 2.2452143271619647e-05     
CG iteration 65, residual = 1.866768382202299e-05     
CG iteration 66, residual = 1.567651790210354e-05     
CG iteration 67, residual = 1.349199170142391e-05     
CG iteration 68, residual = 1.1932267987396287e-05     
CG iteration 69, residual = 1.1006946297625618e-05     
CG iteration 70, residual = 1.0456629339641532e-05     
CG iteration 71, residual = 9.883740543600402e-06     
CG iteration 72, residual = 9.18654477544626e-06     
CG iteration 73, residual = 8.15282914769941e-06     
CG iteration 74, residual = 6.959623023959862e-06     
CG iteration 75, residual = 5.7777984976245495e-06     
CG iteration 76, residual = 4.775007175939122e-06     
CG iteration 77, residual = 3.979591970777657e-06     
CG iteration 78, residual = 3.369765675537563e-06     
CG iteration 79, residual = 2.931468204765155e-06     
CG iteration 80, residual = 2.60634318093901e-06     
CG iteration 81, residual = 2.339109257403547e-06     
CG iteration 82, residual = 2.107326810249669e-06     
CG iteration 83, residual = 1.8783881769549767e-06     
CG iteration 84, residual = 1.6729157949340527e-06     
CG iteration 85, residual = 1.4773593702506955e-06     
CG iteration 86, residual = 1.2787374463305726e-06     
CG iteration 87, residual = 1.1065178086589562e-06     
CG iteration 88, residual = 9.440247565628725e-07     
CG iteration 89, residual = 8.070629830606573e-07     
CG iteration 90, residual = 6.801600792896377e-07     
CG iteration 91, residual = 5.689322461402907e-07     
CG iteration 92, residual = 4.779747968374757e-07     
CG iteration 93, residual = 4.032225331229404e-07     
CG iteration 94, residual = 3.4431548932164006e-07     
CG iteration 95, residual = 2.9836138296028376e-07     
CG iteration 96, residual = 2.648745082625143e-07     
CG iteration 97, residual = 2.383048534478379e-07     
CG iteration 98, residual = 2.131387123002803e-07     
CG iteration 99, residual = 1.898784916121092e-07     
CG iteration 100, residual = 1.6810794227918517e-07     
CG iteration 101, residual = 1.4643428259998588e-07     
CG iteration 102, residual = 1.2650312372658407e-07     
CG iteration 103, residual = 1.0903112531942914e-07     
CG iteration 104, residual = 9.440956877709818e-08     
CG iteration 105, residual = 8.195229678621611e-08     
CG iteration 106, residual = 7.097702472350843e-08     
CG iteration 107, residual = 5.990623454010632e-08     
CG iteration 108, residual = 5.091267024359114e-08     
CG iteration 109, residual = 4.394108542987948e-08     
CG iteration 110, residual = 3.854662919622028e-08     
CG iteration 111, residual = 3.4374907870500334e-08     
CG iteration 112, residual = 3.102602724489258e-08     
CG iteration 113, residual = 2.819916320270096e-08     
CG iteration 114, residual = 2.5958455287589274e-08     
CG iteration 115, residual = 2.377209035003359e-08     
CG iteration 116, residual = 2.154566943384172e-08     
CG iteration 117, residual = 1.8819178101594775e-08     
CG iteration 118, residual = 1.630201406162246e-08     
CG iteration 119, residual = 1.3994410133618132e-08     
CG iteration 120, residual = 1.195825435773463e-08     
CG iteration 121, residual = 1.032613665159515e-08     
CG iteration 122, residual = 9.058935686133658e-09     
CG iteration 123, residual = 8.009682153336291e-09     
CG iteration 124, residual = 7.154863634548925e-09     
CG iteration 125, residual = 6.394862438241238e-09     
CG iteration 126, residual = 5.746973747260621e-09     
CG iteration 127, residual = 5.180094215945222e-09     
CG iteration 128, residual = 4.629842493747649e-09     
CG iteration 129, residual = 4.151238865265757e-09     
CG iteration 130, residual = 3.717168098707984e-09     
CG iteration 131, residual = 3.286366698250553e-09     
CG iteration 132, residual = 2.8759227553543622e-09     
CG iteration 133, residual = 2.4976441033069586e-09     
CG iteration 134, residual = 2.1695613368900353e-09     
CG iteration 135, residual = 1.868479539278659e-09     
CG iteration 136, residual = 1.6186157648527164e-09     
CG iteration 137, residual = 1.3965604961346033e-09     
CG iteration 138, residual = 1.2300628106624675e-09     
CG iteration 139, residual = 1.0837253771781226e-09     
CG iteration 140, residual = 9.728766439066483e-10     
CG iteration 141, residual = 8.732708457756201e-10     
CG iteration 142, residual = 7.863949425366965e-10     
CG iteration 143, residual = 6.947908803722945e-10     
CG iteration 144, residual = 6.07925833981728e-10     
CG iteration 145, residual = 5.24675951257377e-10     
CG iteration 146, residual = 4.50367348582016e-10     
CG iteration 147, residual = 3.8833764485846596e-10     
CG iteration 148, residual = 3.3461607850625935e-10     
CG iteration 149, residual = 2.922510931506465e-10     
CG iteration 150, residual = 2.549684465479315e-10     
CG iteration 151, residual = 2.2459811537180053e-10     
CG iteration 152, residual = 1.9803355549470093e-10     
CG iteration 153, residual = 1.7347919708665725e-10     
CG iteration 154, residual = 1.5224030212052426e-10     
CG iteration 155, residual = 1.335279650796617e-10     
CG iteration 156, residual = 1.1786880569882827e-10     
CG iteration 157, residual = 1.0363448252458666e-10     
CG iteration 158, residual = 9.040220172306008e-11     
CG iteration 159, residual = 7.971559812450795e-11     
CG iteration 160, residual = 7.019193664920986e-11     
CG iteration 161, residual = 6.163245356641879e-11     
CG iteration 162, residual = 5.405299927572636e-11     
CG iteration 163, residual = 4.695454892336523e-11     
CG iteration 164, residual = 4.0594281722006305e-11     
CG iteration 165, residual = 3.4832473581304417e-11     
CG iteration 166, residual = 2.9626569036760813e-11     
CG iteration 167, residual = 2.5116232526283075e-11     
CG iteration 168, residual = 2.1313666734853138e-11     
CG iteration 169, residual = 1.838923944055426e-11     
CG iteration 170, residual = 1.591771140941996e-11     
CG iteration 171, residual = 1.3924783556296984e-11     
CG iteration 172, residual = 1.222737285403413e-11     
CG iteration 173, residual = 1.079455464880582e-11     
CG iteration 174, residual = 9.496967435669054e-12     
CG iteration 175, residual = 8.331575016994831e-12     
CG iteration 176, residual = 7.2199611716040445e-12     
CG iteration 177, residual = 6.2595481293222114e-12     
CG iteration 178, residual = 5.403684066592509e-12     
CG iteration 179, residual = 4.634721546670942e-12     
CG iteration 180, residual = 3.962485904904042e-12     
CG iteration 181, residual = 3.378916537381631e-12     
CG iteration 182, residual = 2.8879516241498604e-12     
CG iteration 183, residual = 2.456948062112392e-12     
CG iteration 184, residual = 2.0874239086725313e-12     
CG iteration 185, residual = 1.793432701881308e-12     
CG iteration 186, residual = 1.5462318450925006e-12     
CG iteration 187, residual = 1.3704534384292498e-12     
CG iteration 188, residual = 1.2300262733780903e-12     
CG iteration 189, residual = 1.1025581842688831e-12     
CG iteration 190, residual = 9.710726594657306e-13     
CG iteration 191, residual = 8.526612361194101e-13     
CG iteration 192, residual = 7.423342942113786e-13     
CG iteration 193, residual = 6.410523312875719e-13     
CG iteration 194, residual = 5.47886358578971e-13     
CG iteration 195, residual = 4.619618411497568e-13     
CG iteration 196, residual = 3.846671392110144e-13     
CG iteration 197, residual = 3.1876724690369686e-13     
CG iteration 198, residual = 2.65851924215303e-13     
CG iteration 199, residual = 2.2442774562309983e-13     
CG iteration 200, residual = 1.9151916963337414e-13     
CG iteration 201, residual = 1.6752304112957501e-13     
CG iteration 202, residual = 1.5074234046264504e-13     
CG iteration 203, residual = 1.3653362423344818e-13     
CG iteration 204, residual = 1.2254065792549437e-13     
CG iteration 205, residual = 1.0945224978282656e-13     
CG iteration 206, residual = 9.857063746295873e-14     
CG iteration 207, residual = 8.823047525152877e-14     
CG iteration 208, residual = 7.812689381991306e-14     
CG iteration 209, residual = 6.678493331640159e-14     
CG iteration 210, residual = 5.5665251475960276e-14     
CG iteration 211, residual = 4.5790061176512394e-14     
CG iteration 212, residual = 3.7889039808345655e-14     
CG iteration 213, residual = 3.167389950326902e-14     
CG iteration 214, residual = 2.6696677388393872e-14     
CG iteration 215, residual = 2.33158312730966e-14     
CG converged in 215 iterations to residual 2.33158312730966e-14
Draw (gfu, draw_vol=False);

Exercise:

  • Experiment with differ problem sizes and preconditioners

  • How big systems can you solve on your computer (watch out for memory usage)